Highly multiplexed imaging acquires the single-cell expression of RNA, metabolites or proteins in a spatially-resolved fashion. These measurements can be visualized across multiple length-scales. First, pixel-level intensities represent the spatial distributions of feature expression with highest resolution. Second, after segmentation, expression values or cell-level metadata (e.g. cell-type information) can be visualized on segmented cell areas. This workflow describes the use of the cytomapper Bioconductor package to demonstrate functions for the visualization of multiplexed read-outs and cell-level information obtained by multiplexed imaging. The main functions of this package allow 1. the visualization of pixel-level information across multiple channels, 2. the display of cell-level information (expression and/or metadata) on segmentation masks and 3. gating and visualisation of single cells.
To follow this tutorial, please visit https://github.com/BodenmillerGroup/cytomapper_demos/tree/main/docs. The compiled .html of this workshop is hosted at: https://bodenmillergroup.github.io/cytomapper_demos. The cytomapper package can be installed via:
if (!requireNamespace("BiocManager", quietly = TRUE))
install.packages("BiocManager")
BiocManager::install("cytomapper")
To reproduce the analysis, clone the repository:
git clone https://github.com/BodenmillerGroup/cytomapper_demos.git
and open the Bioc2021_workshop.Rmd file in the docs folder.
We provide three images, their segmentation masks and the quantified single-cell
data in form of a SingleCellExperiment object in the data folder. The data
is taken from A map of human type 1 diabetes progression by imaging mass
cytometry.
The analysis and visualization of highly-multiplexed imaging data relies on three types of data structures:
Multi-channel images realized as three dimensional arrays or dimensions x,
y, c where the numeric entry of each voxel represents the intensity of the
pixel at position x and y for channel c.
Single-channel segmentation masks realized as matrices where sets of pixels with the same integer ID represent individual objects (here, these are segmented cells)
A table containing the quantified features per cell and channel (e.g. mean pixel intensity)
The cytomapper package handles these data types using objects of the following
S4 classes:
Multiple multi-channel images are stored in form of EBImage::Image objects
within a cytomapper::CytoImageList object
Multiple single-channel segmentation masks are stored in form of
EBImage::Image objects within a cytomapper::CytoImageList object
Per-cell intensity measures and cell-/channel-specific metadata are stored in
a SingleCellExperiment::SingleCellExperiment container
cytomapper overview figure. A) The plotCells function combines a SingleCellExperiment and CytoImageList object to visualize marker expression or cell-specific metadata on segmentation masks. B) The plotPixels function requires a CytoImageList object to visualize the combined expression of up to six markers as composite images
The cytomapper package contains three broad functionalities:
Visualization of pixel-intensities as composites of up to six channels (plotPixels)
Visualization of cell-specific features on segmentation masks (plotCells)
Interactive gating of cells and visualization of gated cells on images (cytomapperShiny)
The follwing demonstration gives an overview on data handling and processing, data visualization and interactive analyses.
We use imaging mass cytometry data
to highlight the functionality of the cytomapper package. However, any imaging
technology is supported as long as the data can be read into R (memory
restrictions and file type restrictions.)
The raw data has been processed using the ImcSegmentationPipeline. In this tutorial, we assume that image segmentation has been performed beforehand.
The following section describes how to read in images (e.g. in .tiff format)
into CytoImageList objects and how to generate SingleCellExperiment
objects from these images.
The cytomapper::loadImages function reads in multi-channel images and segmentation masks into CytoImageList objects.
library(cytomapper)
# Read in 32-bit multi-channel images
(images <- loadImages("../data/images/", pattern = ".tiff"))
## CytoImageList containing 3 image(s)
## names(3): E30_a0_full_clean G23_a0_full_clean J01_a0_full_clean
## Each image contains 38 channel(s)
# Read in 16-bit unsigned integer segmentation masks
(masks <- loadImages("../data/masks/", pattern = ".tiff"))
## CytoImageList containing 3 image(s)
## names(3): E30_a0_full_mask G23_a0_full_mask J01_a0_full_mask
## Each image contains 1 channel
It is always recommended to observe the numeric pixel values to make sure that images were read in correctly:
# multi-channel images - first image, first channel
hist(log10(images[[1]][,,1] + 1))
# Segmentation mask
masks[[1]]
## Image
## colorMode : Grayscale
## storage.mode : double
## dim : 390 452
## frames.total : 1
## frames.render: 1
##
## imageData(object)[1:5,1:6]
## [,1] [,2] [,3] [,4] [,5]
## [1,] 0.000000e+00 0.000000e+00 0.000000e+00 0.000000e+00 4.577707e-05
## [2,] 0.000000e+00 0.000000e+00 0.000000e+00 0.000000e+00 4.577707e-05
## [3,] 4.577707e-05 4.577707e-05 4.577707e-05 4.577707e-05 4.577707e-05
## [4,] 4.577707e-05 4.577707e-05 4.577707e-05 4.577707e-05 4.577707e-05
## [5,] 4.577707e-05 4.577707e-05 4.577707e-05 4.577707e-05 4.577707e-05
## [,6]
## [1,] 0.000000e+00
## [2,] 0.000000e+00
## [3,] 4.577707e-05
## [4,] 4.577707e-05
## [5,] 4.577707e-05
We notice, that the segmentation masks were not read in as integer images.
This behavior arises due to issues with accessing the tiff metadata after
pre-processing using CellProfiler.
In these cases, the cytomapper::scaleImages function can be used to rescale
segmentation masks to only contain integer IDs. Here, the scaling factor
is 2 ^ 16 - 1 = 65535 accounting for the 16-bit unsigned integer encoding.
masks <- scaleImages(masks, 2 ^ 16 - 1)
masks[[1]]
## Image
## colorMode : Grayscale
## storage.mode : double
## dim : 390 452
## frames.total : 1
## frames.render: 1
##
## imageData(object)[1:5,1:6]
## [,1] [,2] [,3] [,4] [,5] [,6]
## [1,] 0 0 0 0 3 0
## [2,] 0 0 0 0 3 0
## [3,] 3 3 3 3 3 3
## [4,] 3 3 3 3 3 3
## [5,] 3 3 3 3 3 3
As an alternative, while reading in the images, the as.is option can be set to TRUE.
masks <- loadImages("../data/masks/", pattern = ".tiff", as.is = TRUE)
masks[[1]]
## Image
## colorMode : Grayscale
## storage.mode : integer
## dim : 390 452
## frames.total : 1
## frames.render: 1
##
## imageData(object)[1:5,1:6]
## [,1] [,2] [,3] [,4] [,5] [,6]
## [1,] 0 0 0 0 3 0
## [2,] 0 0 0 0 3 0
## [3,] 3 3 3 3 3 3
## [4,] 3 3 3 3 3 3
## [5,] 3 3 3 3 3 3
We can already visualize the segmentation mask to get an idea of the tissue structure:
plotCells(masks)
To increase scalability of the cytomapper package, images can be stored on
disk making use of the HDF5Array and DelayedArray package. Reading in the
data can also be parallelised, which is only recommended for large images.
format(object.size(images), units = "Kb")
## [1] "170235.4 Kb"
images_ondisk <- loadImages("../data/images/", pattern = ".tiff",
on_disk = TRUE, h5FilesPath = "../data/images",
BPPARAM = BiocParallel::bpparam())
format(object.size(images_ondisk), units = "Kb")
## [1] "10.1 Kb"
All cytomapper functions support images and masks stored on disk.
SingleCellExperiment objectIn the original analysis of this dataset, single-cell features have been
extracted using CellProfiler. For this tutorial, we provide a
SingleCellExperiment object, which already contains the mean intensities per
cell and channel and all relevant metadata (e.g. cell-type annotation).
The SingleCellExperiment can be read-in from the data folder:
(sce <- readRDS("../data/sce.rds"))
## class: SingleCellExperiment
## dim: 38 5819
## metadata(0):
## assays(2): counts exprs
## rownames(38): H3 SMA ... Ir191 Ir193
## rowData names(15): TubeNb MetalTag ... miCAT2 miCAT
## colnames(5819): E30_1 E30_2 ... J01_2173 J01_2174
## colData names(26): slide id ... Ethnicity BMI
## reducedDimNames(0):
## mainExpName: NULL
## altExpNames(0):
colData(sce)[,c("ImageName", "CellNumber", "CellType")]
## DataFrame with 5819 rows and 3 columns
## ImageName CellNumber CellType
## <character> <integer> <character>
## E30_1 E30 1 ductal
## E30_2 E30 2 macrophage
## E30_3 E30 3 acinar
## E30_4 E30 4 acinar
## E30_5 E30 5 ductal
## ... ... ... ...
## J01_2170 J01 2170 ductal
## J01_2171 J01 2171 acinar
## J01_2172 J01 2172 acinar
## J01_2173 J01 2173 acinar
## J01_2174 J01 2174 acinar
The following section will discuss setting image metadata and image normalization.
Before image visualization, there are a few metadata entries that need to be set.
We will need to set the channel names of the images via the channelNames
getter/setter function. The channel order here is the same as the row order of
the SingleCellExperiment object. We will also need to synchronise the image
IDs across the multi-channel images and segmentation masks by storing a
DataFrame in the elementMetadata slot of the CytoImageList object.
# Add channel names
channelNames(images) <- rownames(sce)
# Add image name to metadata
(mcols(images) <- mcols(masks) <- DataFrame(ImageName = c("E30", "G23", "J01")))
## DataFrame with 3 rows and 1 column
## ImageName
## <character>
## 1 E30
## 2 G23
## 3 J01
Channel normalization is a crucial step for image visualization to enhance
the visibility of biological signals. One option is to transform (e.g. sqrt)
the images. A more widely used alternative is to perform a range-scaling
and clipping of channel intensities.
The cytomapper package exports the normalize function that scales channels
between 0 and 1 across all images (by default). Clipping is performed by setting
the inputRange parameter in the normalize function.
# Min-max normalization
images_norm <- normalize(images)
# Clip the images
images_norm <- normalize(images_norm, inputRange = c(0, 0.2))
The cytomapper::measureObjects function takes the segmentation masks and
multi- channel image objects and computes morphological (e.g. cell shape, size
and location) and intensity features (default: mean intensity per channel and
object/cell).
sce_measured <- measureObjects(mask = masks, image = images,
img_id = "ImageName")
sce_measured
## class: SingleCellExperiment
## dim: 38 5819
## metadata(0):
## assays(1): counts
## rownames(38): H3 SMA ... Ir191 Ir193
## rowData names(0):
## colnames: NULL
## colData names(8): ImageName object_id ... m.majoraxis m.eccentricity
## reducedDimNames(0):
## mainExpName: NULL
## altExpNames(0):
The pixel intensities per cell can be summarized in different ways (e.g. as
quantiles). Furthermore, parallelization is possible by setting BPPARAM = bpparam().
Cell-specific morphological features are stored in colData(sce) while the
mean pixel intensities per cell and channel are stored in counts(sce).
colData(sce_measured)
## DataFrame with 5819 rows and 8 columns
## ImageName object_id s.area s.radius.mean m.cx m.cy
## <character> <numeric> <numeric> <numeric> <numeric> <numeric>
## 1 E30 1 8 1.132373 172.375 1.62500
## 2 E30 2 5 0.969029 150.400 2.60000
## 3 E30 3 85 4.928794 7.400 4.36471
## 4 E30 4 41 3.148247 221.756 3.85366
## 5 E30 5 44 3.454595 190.523 5.95455
## ... ... ... ... ... ... ...
## 5815 J01 2170 40 3.15461 436.950 454.050
## 5816 J01 2171 64 4.26023 457.469 453.688
## 5817 J01 2172 42 3.51294 341.643 454.643
## 5818 J01 2173 10 1.55575 20.900 455.500
## 5819 J01 2174 16 2.46488 44.250 455.625
## m.majoraxis m.eccentricity
## <numeric> <numeric>
## 1 3.96961 0.713215
## 2 4.34170 0.955796
## 3 13.40271 0.767751
## 4 8.17150 0.515845
## 5 10.23176 0.792695
## ... ... ...
## 5815 8.85424 0.722069
## 5816 12.40952 0.812136
## 5817 12.96285 0.936558
## 5818 6.43099 0.928765
## 5819 10.74725 0.985845
counts(sce_measured)[1:10, 1:10]
## [,1] [,2] [,3] [,4] [,5] [,6]
## H3 22.6329811 2.4361224 13.08346010 6.7423978 14.72265212 11.5127982
## SMA 0.7880026 1.0405174 1.26743278 0.9876558 0.32023911 0.7070429
## INS 0.0000000 0.0000000 0.24131852 0.2329416 0.16629238 0.2158641
## CD38 0.3090078 0.2992626 0.32912218 0.1752310 0.20383297 0.2546964
## CD44 2.3301406 0.4281695 0.59169109 0.4615597 1.78500872 0.4622892
## PCSK2 1.1045258 1.0593891 0.56905349 0.5634079 0.85554019 0.7053719
## CD99 0.3115875 0.3936061 0.77623081 1.0573374 0.89659225 0.9886102
## CD68 0.0000000 4.7920089 0.31981109 0.4370047 1.25933963 0.2232831
## MPO 0.0000000 5.1032834 0.07681261 0.1129901 0.04637183 0.1006678
## SLC2A1 2.8200562 2.2636643 0.79776477 1.4814364 1.74246915 0.4607339
## [,7] [,8] [,9] [,10]
## H3 6.65392760 8.7451991 6.6297040 9.66234309
## SMA 1.94122101 1.2825997 1.3894975 0.78360098
## INS 0.09161043 0.1462787 0.1709020 0.09961374
## CD38 0.26884190 0.5762133 0.3820182 0.50709601
## CD44 0.40098510 0.5384129 0.3380702 0.75498592
## PCSK2 0.65091868 0.8107124 1.0038303 0.68700536
## CD99 0.73302101 0.7973502 0.5520242 0.94204760
## CD68 1.04103449 0.7391340 0.2099665 1.36866162
## MPO 3.35050052 5.6631805 0.6724382 8.20452772
## SLC2A1 1.61920737 2.3110813 1.0785549 2.68912075
After having formatted the data, we can now move on to data visualization.
The cytomapper::plotPixels function is the main function for visualization
of pixel intensities.
The easiest way of visualizing images is to just specify channel names:
plotPixels(images,
colour_by = c("PIN", "CDH"))
Here, the intensities appear weak. For this, we have normalized the images:
plotPixels(images_norm,
colour_by = c("PIN", "CDH"))
We can see a decline in pro-insulin expression while E-cadherin expression stays the same.
Instead of normalizing the images, you can also specify the background (b), contrast (c) and gamma (g) levels of the images:
plotPixels(images,
colour_by = c("PIN", "CDH"),
bcg = list(PIN = c(0, 5, 1),
CDH = c(0, 5, 1)))
The cytomapper package allows flexible adjustments of all visual attributes of the images:
plotPixels(
image = images,
colour_by = c("PIN", "CD4", "CD8a"),
colour = list(PIN = c("black", "yellow"),
CD4 = c("black", "blue"),
CD8a = c("black", "red")),
bcg = list(PIN = c(0, 10, 1),
CD4 = c(0, 8, 1),
CD8a = c(0, 10, 1)),
image_title = list(
text = c("Non-diabetic",
"Recent onset T1D",
"Long duration T1D")
),
scale_bar = list(
length = 100,
label = expression("100 " ~ mu * "m")
))
The colour_by parameter defines the channel names by which to colour the
composite. Per channel, a colour scale is generated by setting colour. The
attributes of the image titles can be set via the parameter image_title and
attributes of the scale bar are set via the parameter scale_bar.
To see all available parameter options, refer to ?plotting-param
In the next step, we can visualize mean pixel intensities and cell-specific metadata directly on the segmentation masks.
For this, we will use the cytomapper::plotCells function while combining the
SingleCellExperiment object and the CytoImageList storing the segmentation
masks:
plotCells(mask = masks, object = sce,
cell_id = "CellNumber", img_id = "ImageName",
colour_by = c("PIN", "CD8a", "CD4"),
exprs_values = "exprs")
The assay(sce, "exprs") slot stores the asinh-transformed mean pixel intensities
per channel and cell. The cell_id entry needs to be stored in the colData(sce)
slot to link each individual cell to their corresponding pixels in the segmentation mask.
The img_id parameter specifies the colData(sce) entry storing the name/ID
of the image to which each cell belongs. This should also be an entry in
mcols(masks).
The SingleCellExperiment object is subsettable to highlight only certain cells:
cur_sce <- sce[,sce$CellType %in%
c("beta", "alpha", "delta", "Tc", "Th")]
plotCells(mask = masks, object = cur_sce,
cell_id = "CellNumber", img_id = "ImageName",
colour_by = c("PIN", "CD8a", "CD4"),
exprs_values = "exprs",
missing_colour = "white")
Here, missing_colour defines the colour of cells missing in the SingleCellExperiment object.
This can also be useful when visualizing cell-specific metadata:
plotCells(
mask = masks,
object = cur_sce,
cell_id = "CellNumber",
img_id = "ImageName",
colour_by = "CellType",
image_title = list(
text = c("Non-diabetic",
"Recent onset T1D",
"Long duration T1D"),
colour = "black"),
scale_bar = list(
length = 100,
label = expression("100 " ~ mu * "m"),
colour = "black"),
missing_colour = "white",
background_colour = "gray")
Finally, all three objects (multi-channel images, segmentation masks and the single-cell data) can be combined to outline cells on composite images.
As a first demonstration, we select only pancreatic beta cells and visualize their corresponing marker pro-insulin:
cur_sce <- sce[,sce$CellType == "beta"]
plotPixels(images_norm,
mask = masks,
object = cur_sce,
cell_id = "CellNumber",
img_id = "ImageName",
colour_by = c("PIN", "H3"),
bcg = list(H3 = c(0, 2, 1)),
colour = list(PIN = c("black", "yellow"),
H3 = c("black", "blue"),
CellType = c(beta = "red")),
outline_by = "CellType")
One can now also visualize multiple cell types and multiple markers together:
We have developed a shiny application to gate cells based on their mean
pixel intensity with additional visualization on the images.
The idea behind this gating strategy is to generate ground truth cell type labels, which can be used to train a classifier and classify cells rather than relying on clustering.
The shiny app can be opened using the cytomapper::cytomapperShiny function:
if (interactive()) {
cytomapperShiny(sce, mask = masks, image = images,
img_id = "ImageName", cell_id = "CellNumber")
}
Using the Shiny application
The help page provides a recommended workflow on how to most efficiently use the
app. The workflow is solely a recommendation - the app provides full flexibility
to change settings during each step. To see the full documentation, please refer
to the help page found at ?cytomapperShiny
Select the number of plots The slider under “General controls” can be used to specify the number of plots on which to perform gating. Up to two markers can be visualized per plot.
Select the sample
The assay dropdown selection under “General controls” allows the user to specify
on which assay entry to perform gating. In most cases, a log- or
arcsinh-transformation can help to distinguish between ‘positive’ and ‘negative’
populations.
Select the markers
For each plot, up to two markers can be specified. If selecting a single marker,
please specify this marker in the first of the two dropdown menus. A violin plot
is used to visualize the expression of a single marker while a scatter plot is
used to visualize the expression of two markers.
Gate cells
When selecting cells in one plot, only those cells are visualized on the
following plot. Once markers, the assay or the number of plots are changed,
gates are cleared.
Observe the selected cells
After gating, the selected cells are visualized on the corresponding images by
switching to the “Images tab”. By default, the first marker is selected. The
user can change the displayed marker or press reset marker to switch to the
markers used for gating. If a multi-channel image object is provided, the
contrast of the image can be changed. The right panel visualizes the selected
cells either by filling in the segmentation masks or by outlining the cells on
the images.
Change samples
Samples can now be iteratively changed using the dropdown menu under General
controls . The gates will remain on the plots and can be adjusted for each
sample.
Save the selected cells
Finally, the selected cells can be saved by clicking the download button next to
the ‘?’ symbol. The selected cells will be stored as a SingleCellExperiment
object in .rds format. Per selection, the user can provide a Cell label that
will be stored in the colData under the cytomapper_CellLabel entry of the
downloaded object.
For pre-processing multiplexed imaging data, please refer to the ImcSegmentationPipeline and/or the steinbock package.
To test the cytomapper package on different datasets, check out the imcdatasets Bioconductor package.
The imcRtools package is currently being developed to facilitate handling of multiplexed imaging data and spatial analysis.
## R version 4.1.0 (2021-05-18)
## Platform: x86_64-apple-darwin17.0 (64-bit)
## Running under: macOS Catalina 10.15.7
##
## Matrix products: default
## BLAS: /Library/Frameworks/R.framework/Versions/4.1/Resources/lib/libRblas.dylib
## LAPACK: /Library/Frameworks/R.framework/Versions/4.1/Resources/lib/libRlapack.dylib
##
## locale:
## [1] en_US.UTF-8/en_US.UTF-8/en_US.UTF-8/C/en_US.UTF-8/en_US.UTF-8
##
## attached base packages:
## [1] parallel stats4 stats graphics grDevices utils datasets
## [8] methods base
##
## other attached packages:
## [1] cytomapper_1.4.1 SingleCellExperiment_1.14.1
## [3] SummarizedExperiment_1.22.0 Biobase_2.52.0
## [5] GenomicRanges_1.44.0 GenomeInfoDb_1.28.1
## [7] IRanges_2.26.0 S4Vectors_0.30.0
## [9] BiocGenerics_0.38.0 MatrixGenerics_1.4.0
## [11] matrixStats_0.60.0 EBImage_4.34.0
## [13] BiocStyle_2.20.2
##
## loaded via a namespace (and not attached):
## [1] bitops_1.0-7 fontawesome_0.2.2 RColorBrewer_1.1-2
## [4] tools_4.1.0 bslib_0.2.5.1 utf8_1.2.2
## [7] R6_2.5.0 svgPanZoom_0.3.4 HDF5Array_1.20.0
## [10] vipor_0.4.5 DBI_1.1.1 colorspace_2.0-2
## [13] rhdf5filters_1.4.0 raster_3.4-13 sp_1.4-5
## [16] gridExtra_2.3 tidyselect_1.1.1 compiler_4.1.0
## [19] DelayedArray_0.18.0 bookdown_0.22 sass_0.4.0
## [22] scales_1.1.1 systemfonts_1.0.2 stringr_1.4.0
## [25] digest_0.6.27 tiff_0.1-8 fftwtools_0.9-11
## [28] rmarkdown_2.9 svglite_2.0.0 XVector_0.32.0
## [31] jpeg_0.1-9 pkgconfig_2.0.3 htmltools_0.5.1.1
## [34] highr_0.9 fastmap_1.1.0 htmlwidgets_1.5.3
## [37] rlang_0.4.11 shiny_1.6.0 jquerylib_0.1.4
## [40] generics_0.1.0 jsonlite_1.7.2 BiocParallel_1.26.1
## [43] dplyr_1.0.7 RCurl_1.98-1.3 magrittr_2.0.1
## [46] GenomeInfoDbData_1.2.6 Matrix_1.3-3 Rcpp_1.0.7
## [49] ggbeeswarm_0.6.0 munsell_0.5.0 Rhdf5lib_1.14.2
## [52] fansi_0.5.0 viridis_0.6.1 abind_1.4-5
## [55] lifecycle_1.0.0 stringi_1.7.3 yaml_2.2.1
## [58] zlibbioc_1.38.0 rhdf5_2.36.0 grid_4.1.0
## [61] promises_1.2.0.1 shinydashboard_0.7.1 crayon_1.4.1
## [64] lattice_0.20-44 magick_2.7.2 locfit_1.5-9.4
## [67] knitr_1.33 pillar_1.6.1 codetools_0.2-18
## [70] glue_1.4.2 evaluate_0.14 BiocManager_1.30.16
## [73] png_0.1-7 vctrs_0.3.8 httpuv_1.6.1
## [76] gtable_0.3.0 purrr_0.3.4 assertthat_0.2.1
## [79] ggplot2_3.3.5 xfun_0.24 mime_0.11
## [82] xtable_1.8-4 later_1.2.0 viridisLite_0.4.0
## [85] tibble_3.1.3 beeswarm_0.4.0 ellipsis_0.3.2